home *** CD-ROM | disk | FTP | other *** search
/ Mac Magazin/MacEasy 21 / Mac Magazin and MacEasy Magazine CD - Issue 21.iso / Wissenschaft & Technik / yorick12vr1-nofpu folder / startup / matrix.i < prev    next >
Text File  |  1995-07-26  |  17KB  |  497 lines

  1. /*
  2.    MATRIX.I
  3.    Yorick interface to LAPACK matrix solving routines.
  4.  
  5.    $Id: matrix.i,v 1.1 1993/08/07 19:57:08 munro Exp $
  6.  */
  7. /*    Copyright (c) 1994.  The Regents of the University of California.
  8.                     All rights reserved.  */
  9.  
  10. /* ------------------------------------------------------------------------ */
  11.  
  12. func unit(n, m)
  13. /* DOCUMENT unit(n)
  14.          or unit(n, m)
  15.      returns n-by-n (or n-by-m) unit matrix, i.e.- matrix with diagonal
  16.      elements all 1.0, off diagonal elements 0.0
  17.  */
  18. {
  19.   if (is_void(m)) m= n;
  20.   u= array(0.0, n, m);
  21.   u(1:numberof(u):n+1)= 1.0;
  22.   return u;
  23. }
  24.  
  25. /* ------------------------------------------------------------------------ */
  26.  
  27. func TDsolve(c, d, e, b, which=)
  28. /* DOCUMENT TDsolve(c, d, e, b)
  29.          or TDsolve(c, d, e, b, which=which)
  30.  
  31.      returns the solution to the tridiagonal system:
  32.         D(1)*x(1)       + E(1)*x(2)                       = B(1)
  33.     C(1:-1)*x(1:-2) + D(2:-1)*x(2:-1) + E(2:0)*x(3:0) = B(2:-1)
  34.                       C(0)*x(-1)      + D(0)*x(0)     = B(0)
  35.      (i.e.- C is the subdiagonal, D the diagonal, and E the superdiagonal;
  36.      C and E have one fewer element than D, which is the same length as
  37.      both B and x)
  38.  
  39.      B may have additional dimensions, in which case the returned x
  40.      will have the same additional dimensions.  The WHICH dimension of B,
  41.      and of the returned x is the one of length n which participates
  42.      in the matrix solve.  By default, WHICH=1, so that the equations
  43.      being solved involve B(,..) and x(+,..).
  44.      Non-positive WHICH counts from the final dimension (as for the
  45.      sort and transpose functions), so that WHICH=0 involves B(..,)
  46.      and x(..,+).
  47.  
  48.      The C, D, and E arguments may be either scalars or vectors; they
  49.      will be broadcast as appropriate.
  50.  
  51.   SEE ALSO: LUsolve, QRsolve, SVsolve, SVdec
  52.  */
  53. {
  54.   /* check validity of b argument */
  55.   if (structof(b)==complex) error, "expecting a non-complex RHS vector";
  56.   dims= dimsof(b);
  57.   ndb= is_void(dims)? 0 : dims(1);
  58.   if (is_void(which)) which= 1;
  59.   else if (which<=0) which+= ndb;
  60.   if (!ndb) error, "RHS must have at least one dimension";
  61.   n= dims(1+which);
  62.   b= double(b);   /* copy of RHS to be transformed into solution */
  63.   nrhs= numberof(b)/n;
  64.  
  65.   /* put first matrix dimension of b first */
  66.   if (which!=1) b= transpose(b, [1,which]);
  67.  
  68.   /* copy, force to double, and broadcast matrix diagonals
  69.      -- also will blow up on conformability error */
  70.   cc= ee= array(0.0, n-1);
  71.   dd= array(0.0, n);
  72.   cc()= c;
  73.   dd()= d;
  74.   ee()= e;
  75.  
  76.   info= 0;
  77.   _dgtsv, n, nrhs, cc, dd, ee, b, n, info;
  78.   if (info) error, "tridiagonal element "+pr1(info)+" of became 0.0";
  79.  
  80.   /* restore proper order of result if necessary */
  81.   if (which!=1) b= transpose(b, [1,which]);
  82.  
  83.   return b;
  84. }
  85.  
  86. extern _dgtsv;
  87. /* PROTOTYPE
  88.    void dgtsv(long n, long nrhs, double array c, double array d,
  89.               double array e, double array b, long ldb, long array info)
  90.  */
  91. /* DOCUMENT _dgtsv
  92.      LAPACK dgtsv routine.
  93.  */
  94.  
  95. /* ------------------------------------------------------------------------ */
  96.  
  97. func LUsolve(a, b, which=, job=, det=)
  98. /* DOCUMENT LUsolve(a, b)
  99.          or LUsolve(a, b, which=which)
  100.      or a_inverse= LUsolve(a)
  101.  
  102.      returns the solution x of the matrix equation:
  103.         A(,+)*x(+) = B
  104.      If A is an n-by-n matrix then B must have length n, and the returned
  105.      x will also have length n.
  106.  
  107.      B may have additional dimensions, in which case the returned x
  108.      will have the same additional dimensions.  The WHICH dimension of B,
  109.      and of the returned x is the one of length n which participates
  110.      in the matrix solve.  By default, WHICH=1, so that the equations
  111.      being solved are:
  112.         A(,+)*x(+,..) = B
  113.      Non-positive WHICH counts from the final dimension (as for the
  114.      sort and transpose functions), so that WHICH=0 solves:
  115.         x(..,+)*A(,+) = B
  116.  
  117.      If the B argument is omitted, the inverse of A is returned:
  118.      A(,+)*x(+,) and A(,+)*x(,+) will be unit matrices.
  119.  
  120.      LUsolve works by LU decomposition using Gaussian elimination with
  121.      pivoting.  It is the fastest way to solve square matrices.  QRsolve
  122.      handles non-square matrices, as does SVsolve.  SVsolve is slowest,
  123.      but can deal with highly singular matrices sensibly.
  124.  
  125.    SEE ALSO: QRsolve, TDsolve, SVsolve, SVdec, LUrcond
  126.  */
  127. {
  128.   /* get n, m, dims, nrhs, checking validity of a and b */
  129.   { local dims, n, m, nrhs; }
  130.   _get_matrix, 1;
  131.   if (m!=n) error, "expecting a square matrix";
  132.  
  133.   if (is_void(b)) {
  134.     b= unit(n);
  135.     nrhs= n;
  136.     which= 1;
  137.   }
  138.  
  139.   /* perform LU solve */
  140.   pivot= array(0, n);
  141.   info= 0;
  142.   _dgesv, n, nrhs, a, n, pivot, b, n, info;
  143.   /* row i interchanged with row pivot(i) --> permutation matrix P
  144.      a now contains the L and U factors from the decomposition;
  145.      original a= P*L*U */
  146.   if (info) error, "matrix is (numerically) singular";
  147.  
  148.   /* restore proper order of result if necessary */
  149.   if (which!=1) b= transpose(b, [1,which]);
  150.  
  151.   return b;
  152. }
  153.  
  154. func LUrcond(a, one_norm=)
  155. /* DOCUMENT LUrcond(a)
  156.          or LUrcond(a, one_norm=1)
  157.      returns the reciprocal condition number of the N-by-N matrix A.
  158.      If the ONE_NORM argument is non-nil and non-zero, the 1-norm
  159.      condition number is returned, otherwise the infinity-norm condition
  160.      number is returned.
  161.  
  162.      The condition number is the ratio of the largest to the smallest
  163.      singular value, max(singular_values)*max(1/singular_values) (or
  164.      sum(abs(singular_values)*sum(abs(1/singular_values)) if ONE_NORM
  165.      is selected?).  If the reciprocal condition number is near zero
  166.      then A is numerically singular; specifically, if
  167.           1.0+LUrcond(a) == 1.0
  168.      then A is numerically singular.
  169.  
  170.    SEE ALSO: LUsolve
  171.  */
  172. {
  173.   dims= dimsof(a);
  174.   if (is_void(dims) || dims(1)!=2 || dims(2)!=dims(3) ||
  175.       structof(a)==complex)
  176.     error, "expecting a square 2D real matrix";
  177.   n= dims(2);
  178.   a= double(a);
  179.   pivot= array(0, n);
  180.   info= 0;
  181.   _dgetrf, n, n, a, n, pivot, info;
  182.   work= array(double, 4*n);
  183.   iwork= array(0, n);
  184.   rcond= 0.0;
  185.   if (!one_norm) {
  186.     one_norm= 0;
  187.     anorm= abs(a)(max,sum);
  188.   } else {
  189.     one_norm= 1;
  190.     anorm= abs(a)(sum,max);
  191.   }
  192.   _dgecox, one_norm, n, a, n, anorm, rcond, work, iwork, info;
  193.   return rcond;
  194. }
  195.  
  196. extern _dgesv;
  197. /* PROTOTYPE
  198.    void dgesv(long n, long nrhs, double array a, long lda,
  199.               long array pivot, double array b, long ldb, long array info)
  200.  */
  201. /* DOCUMENT _dgesv
  202.      LAPACK dgesv routine.
  203.  */
  204.  
  205. extern _dgetrf;
  206. /* PROTOTYPE
  207.    void dgetrf(long m, long n, double array a, long lda,
  208.                long array pivot, long array info)
  209.  */
  210. /* DOCUMENT _dgetrf
  211.      LAPACK dgetrf routine.  Performs LU factorization.
  212.  */
  213.  
  214. extern _dgecox;
  215. /* PROTOTYPE
  216.    void dgecox(long norm, long n, double array a, long lda,
  217.                double anorm, double array rcond, double array work,
  218.                long array iwork, long array info)
  219.  */
  220. /* DOCUMENT _dgecox
  221.      LAPACK dgecon routine, except norm argument not a string.
  222.  */
  223.  
  224. /* ------------------------------------------------------------------------ */
  225.  
  226. func QRsolve(a, b, which=)
  227. /* DOCUMENT QRsolve(a, b)
  228.          or QRsolve(a, b, which=which)
  229.  
  230.      returns the solution x (in a least squares or least norm sense
  231.      described below) of the matrix equation:
  232.         A(,+)*x(+) = B
  233.      If A is an m-by-n matrix (i.e.- m equations in n unknowns), then B
  234.      must have length m, and the returned x will have length n.
  235.  
  236.      If n<m, the system is overdetermined -- no solutions are possible
  237.              -- the returned x minimizes sqrt(sum((A(,+)*x(+) - B)^2))
  238.      If n>m, the system is underdetermined -- many solutions are possible
  239.              -- the returned x has minimum L2 norm among all solutions
  240.  
  241.      B may have additional dimensions, in which case the returned x
  242.      will have the same additional dimensions also have those dimensions.
  243.      The WHICH dimension of B and the returned x is the one of length m
  244.      or n which participates in the matrix solve.  By default, WHICH=1,
  245.      so that the equations being solved are:
  246.         A(,+)*x(+,..) = B
  247.      Non-positive WHICH counts from the final dimension (as for the
  248.      sort and transpose functions), so that WHICH=0 solves:
  249.         A(,+)*x(..,+) = B
  250.  
  251.      QRsolve works by QR factorization if n<m, or LQ factorization if n>m.
  252.      QRsolve is slower than LUsolve.  Its main attraction is that it can
  253.      handle overdetermined or underdetermined systems of equations
  254.      (nonsquare matrices).  QRsolve may fail for singular systems; try
  255.      SVsolve in this case.
  256.  
  257.    SEE ALSO: LUsolve, TDsolve, SVsolve, SVdec
  258.  */
  259. {
  260.   /* get n, m, dims, nrhs, checking validity of a and b */
  261.   { local dims, n, m, nrhs; }
  262.   _get_matrix, 0;
  263.  
  264.   /* set up and perform QR or LQ solve --
  265.      first call returns optimal workspace length */
  266.   work= 0.0;
  267.   info= 0;
  268.   mnmax= max(m,n);
  269.   _dgelx, 0, m, n, nrhs, a, m, b, mnmax, work, 1, info;
  270.   if (info==-10) {
  271.     lwork= long(work);
  272.     work= array(0.0, lwork);
  273.     _dgelx, 0, m, n, nrhs, a, m, b, mnmax, work, lwork, info;
  274.   }
  275.   if (info) error, "matrix is (numerically) singular"; /* impossible? */
  276.  
  277.   /* restore proper order of result if necessary */
  278.   if (n<mnmax) b= b(1:n,..)
  279.   if (which!=1) b= transpose(b, [1,which]);
  280.  
  281.   return b;
  282. }
  283.  
  284. extern _dgelx;
  285. /* PROTOTYPE
  286.    void dgelx(long trans, long m, long n, long nrhs,
  287.               double array a, long lda, double array b, long ldb,
  288.               double array work, long lwork, long array info)
  289.  */
  290. /* DOCUMENT _dgelx
  291.      LAPACK dgels routine, except trans argument not a string.
  292.  */
  293.  
  294. /* ------------------------------------------------------------------------ */
  295.  
  296. func SVsolve(a, b, rcond, which=)
  297. /* DOCUMENT SVsolve(a, b)
  298.          or SVsolve(a, b, rcond)
  299.          or SVsolve(a, b, rcond, which=which)
  300.  
  301.      returns the solution x (in a least squares sense described below) of
  302.      the matrix equation:
  303.         A(,+)*x(+) = B
  304.      If A is an m-by-n matrix (i.e.- m equations in n unknowns), then B
  305.      must have length m, and the returned x will have length n.
  306.  
  307.      If n<m, the system is overdetermined -- no solutions are possible
  308.              -- the returned x minimizes sqrt(sum((A(,+)*x(+) - B)^2))
  309.      If n>m, the system is underdetermined -- many solutions are possible
  310.              -- the returned x has minimum L2 norm among all solutions
  311.  
  312.      SVsolve works by singular value decomposition, therefore it is
  313.      immune to failure due to singularity of the A matrix.  The optional
  314.      RCOND argument defaults to 1.0e-9; singular values less than RCOND
  315.      times the largest singular value (absolute value) will be set to 0.0.
  316.      If RCOND<=0.0, machine precision is used.  The effective rank of the
  317.      matrix is returned as the external variable SVrank.
  318.  
  319.      You can examine the details of the SVD by calling the SVdec routine,
  320.      which returns the singular vectors as well as the singular values.
  321.      Numerical Recipes (Press, et. al. Cambridge University Press 1988)
  322.      has a good discussion of how to use the SVD -- see section 2.9.
  323.  
  324.      B may have additional dimensions, in which case the returned x
  325.      will have the same additional dimensions.  The WHICH argument
  326.      (default 1) controls which dimension of B takes part in the matrix
  327.      solve.  See QRsolve or LUsolve for a complete discussion.
  328.  
  329.    SEE ALSO: SVdec, LUsolve, QRsolve, TDsolve
  330.  */
  331. {
  332.   /* get n, m, dims, nrhs, checking validity of a and b */
  333.   { local dims, n, m, nrhs; }
  334.   _get_matrix, 0;
  335.  
  336.   if (is_void(rcond)) rcond= 1.e-9;
  337.   else rcond= double(rcond);
  338.  
  339.   /* set up and perform SVD solve --
  340.      first call returns optimal workspace length */
  341.   { extern SVrank; }
  342.   work= 0.0;
  343.   info= SVrank= 0;
  344.   s= array(0.0, min(m,n));
  345.   mnmax= max(m,n);
  346.   _dgelss, m, n, nrhs, a, m, b, mnmax, s, rcond, SVrank, work, 1, info;
  347.   if (info==-12) {
  348.     lwork= long(work);
  349.     work= array(0.0, lwork);
  350.     _dgelss, m, n, nrhs, a, m, b, mnmax, s, rcond, SVrank, work, lwork, info;
  351.   }
  352.   if (info) error, "SVD algorithm failed to converge - wow";
  353.  
  354.   /* restore proper order of result if necessary */
  355.   if (n<mnmax) b= b(1:n,..)
  356.   if (which!=1) b= transpose(b, [1,which]);
  357.  
  358.   return b;
  359. }
  360.  
  361. func SVdec(a, &u, &vt, full=)
  362. /* DOCUMENT s= SVdec(a, u, vt)
  363.          or s= SVdec(a, u, vt, full=1)
  364.  
  365.      performs the singular value decomposition of the m-by-n matrix A:
  366.         A = (U(,+) * SIGMA(+,))(,+) * VT(,+)
  367.      where U is an m-by-m orthogonal matrix, VT is an n-by-n orthogonal
  368.      matrix, and SIGMA is an m-by-n matrix which is zero except for its
  369.      min(m,n) diagonal elements.  These diagonal elements are the return
  370.      value of the function, S.  The returned S is always arranged in
  371.      order of descending absolute value.  U(,1:min(m,n)) are the left
  372.      singular vectors corresponding to the min(m,n) elements of S;
  373.      VT(1:min(m,n),) are the right singular vectors.  (The original A
  374.      matrix maps a right singular vector onto the corresponding left
  375.      singular vector, stretched by a factor of the singular value.)
  376.  
  377.      Note that U and VT are strictly outputs; if you don't need them,
  378.      they need not be present in the calling sequence.
  379.  
  380.      By default, U will be an m-by-min(m,n) matrix, and V will be
  381.      a min(m,n)-by-n matrix (i.e.- only the singular vextors are returned,
  382.      not the full orthogonal matrices).  Set the FULL keyword to a
  383.      non-zero value to get the full m-by-m and n-by-n matrices.
  384.  
  385.      On rare occasions, the routine may fail; if it does, the
  386.      first SVinfo values of the returned S are incorrect.  Hence,
  387.      the external variable SVinfo will be 0 after a successful call
  388.      to SVdec.  If SVinfo>0, then external SVe contains the superdiagonal
  389.      elements of the bidiagonal matrix whose diagonal is the returned
  390.      S, and that bidiagonal matrix is equal to (U(+,)*A(+,))(,+) * V(+,).
  391.  
  392.      Numerical Recipes (Press, et. al. Cambridge University Press 1988)
  393.      has a good discussion of how to use the SVD -- see section 2.9.
  394.  
  395.    SEE ALSO: SVsolve, LUsolve, QRsolve, TDsolve
  396.  */
  397. {
  398.   /* get n, m, dims, nrhs, checking validity of a and b */
  399.   { local dims, n, m, nrhs; }
  400.   b= [];
  401.   _get_matrix, 1;
  402.  
  403.   if (!full) full= 0;
  404.   else full= 1;
  405.  
  406.   /* set up and perform SVD solve --
  407.      first call returns optimal workspace length */
  408.   work= 0.0;
  409.   info= 0;
  410.   s= array(0.0, min(m,n));
  411.   if (full) {
  412.     u= array(0.0, m, m);
  413.     vt= array(0.0, n, n);
  414.     ldvt= n;
  415.   } else {
  416.     ldvt= min(m, n);
  417.     u= array(0.0, m, ldvt);
  418.     vt= array(0.0, ldvt, n);
  419.   }
  420.   _dgesvx, full, m, n, a, m, s, u, m, vt, ldvt, work, 1, info;
  421.   if (info==-13) {
  422.     lwork= long(work);
  423.     work= array(0.0, lwork);
  424.     _dgesvx, full, m, n, a, m, s, u, m, vt, ldvt, work, lwork, info;
  425.   }
  426.   if (info) error, "SVD algorithm failed to converge - wow";
  427.  
  428.   return s;
  429. }
  430.  
  431. extern _dgelss;
  432. /* PROTOTYPE
  433.    void dgelss(long m, long n, long nrhs, double array a, long lda,
  434.                double array b, long ldb, double array s, double rcond,
  435.            long array rank, double array work, long lwork,
  436.            long array info)
  437.  */
  438. /* DOCUMENT _dgelss
  439.      LAPACK dgelss routine.
  440.  */
  441.  
  442. extern _dgesvx;
  443. /* PROTOTYPE
  444.    void dgesvx(long job, long m, long n, double array a, long lda,
  445.                double array s, double array u, long ldu,
  446.                double array vt, long ldvt, double array work, long lwork,
  447.            long array info)
  448.  */
  449. /* DOCUMENT _dgesvx
  450.      LAPACK dgesvd routine, except jobu and jobvt are not strings.
  451.  */
  452.  
  453. /* ------------------------------------------------------------------------ */
  454.  
  455. func _get_matrix(b_optional)
  456. {
  457.   { extern dims, n, m, nrhs; }
  458.  
  459.   /* check validity of a argument */
  460.   dims= dimsof(a);
  461.   if (dims(1)!=2 || structof(a)==complex)
  462.     error, "expecting a non-complex 2D matrix";
  463.   a= double(a);  /* copy a to avoid clobbering, as well as force type */
  464.   m= dims(2);
  465.   n= dims(3);
  466.  
  467.   /* check validity of b argument */
  468.   if (!b_optional || !is_void(b)) {
  469.     dims= dimsof(b);
  470.     ndb= is_void(dims)? 0 : dims(1);
  471.     if (is_void(which)) which= 1;
  472.     else if (which<=0) which+= ndb;
  473.     if (!ndb || dims(1+which)!=m)
  474.       error, "RHS dimensions not conformable with matrix dimensions";
  475.     if (structof(b)==complex) error, "expecting a non-complex RHS vector";
  476.     b= double(b);  /* copy to avoid clobbering, and force type */
  477.     nrhs= numberof(b)/m;
  478.  
  479.     /* put first matrix dimension of b first */
  480.     if (which!=1) b= transpose(b, [1,which]);
  481.  
  482.     /* be sure that the first dimension of b is at least n */
  483.     if (n>m) {
  484.       dims= dimsof(b);
  485.       dims(2)= n;
  486.       bn= array(0.0, dims);
  487.       bn(1:m,..)= b;
  488.       b= bn;
  489.     }
  490.  
  491.   } else {
  492.     nrhs= 0;
  493.   }
  494. }
  495.  
  496. /* ------------------------------------------------------------------------ */
  497.